# Copyright (c) HySoP 2011-2024
#
# This file is part of HySoP software.
# See "https://particle_methods.gricad-pages.univ-grenoble-alpes.fr/hysop-doc/"
# for further info.
#
# Licensed under the Apache License, Version 2.0 (the "License");
# you may not use this file except in compliance with the License.
# You may obtain a copy of the License at
#
# http://www.apache.org/licenses/LICENSE-2.0
#
# Unless required by applicable law or agreed to in writing, software
# distributed under the License is distributed on an "AS IS" BASIS,
# WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
# See the License for the specific language governing permissions and
# limitations under the License.
"""
@file min_max.py
MinMaxFieldStatisticsBase: compute min(F), max(F) and/or max(|F|) for a given field
MinMaxGradientStatisticsBase: compute min(gradF), max(gradF) and/or max(|gradF|) for a given field, direction-wise.
"""
import numpy as np
from abc import abstractmethod
from hysop.tools.htypes import check_instance, first_not_None, to_tuple
from hysop.tools.henum import EnumFactory
from hysop.tools.decorators import debug
from hysop.fields.continuous_field import Field
from hysop.constants import Backend
from hysop.parameters.tensor_parameter import TensorParameter
from hysop.topology.cartesian_descriptor import CartesianTopologyDescriptors
from hysop.core.graph.computational_operator import ComputationalGraphOperator
from hysop.parameters.scalar_parameter import ScalarParameter
from hysop.core.mpi import MPI
[docs]
class MinMaxFieldStatisticsBase:
"""
Abstract operator base to compute min and max statistics on a specific field.
"""
[docs]
@classmethod
def supports_multiple_topologies(cls):
return True
[docs]
@classmethod
def build_parameters(
cls,
field,
components,
all_quiet,
Fmin,
Fmax,
Finf,
pbasename,
ppbasename,
dtype=None,
):
if (
((Fmin is None) or (Fmin is False))
and ((Fmax is None) or (Fmax is False))
and ((Finf is None) or (Finf is False))
):
msg = "No statistics were requested."
msg += "\nPlease specify Fmin, Fmax and/or Finf by either setting "
msg += " their value to True, or by by passing an already existing "
msg += " tensor parameter."
raise ValueError(msg)
if field is not None:
dtype = first_not_None(dtype, field.dtype)
pbasename = first_not_None(pbasename, field.name)
ppbasename = first_not_None(ppbasename, field.pretty_name)
def make_param(k, quiet):
return TensorParameter(
name=names[k],
pretty_name=pretty_names[k],
dtype=field.dtype,
shape=(nb_components,),
quiet=quiet,
)
names = {
"Fmin": f"{pbasename}min",
"Fmax": f"{pbasename}max",
"Finf": f"{pbasename}inf",
}
pretty_names = {
"Fmin": f"{ppbasename}₋",
"Fmax": f"{ppbasename}₊",
"Finf": f"|{ppbasename}|₊",
}
if field is not None:
components = first_not_None(components, range(field.nb_components))
components = to_tuple(components)
nb_components = len(components)
parameters = {}
_parameters = dict(Fmin=Fmin, Fmax=Fmax, Finf=Finf)
for k, v in _parameters.items():
param = _parameters[k]
if isinstance(param, TensorParameter):
pass
elif param is True:
param = make_param(k, quiet=all_quiet)
elif (_parameters["Finf"] is not None) and ((k == "Fmin") or (k == "Fmax")):
param = make_param(k, quiet=True)
else:
param = None
if param is not None:
assert np.prod(param.shape) == nb_components
parameters[k] = param
return parameters
@debug
def __new__(
cls,
field,
components=None,
coeffs=None,
Fmin=None,
Fmax=None,
Finf=None,
all_quiet=None,
name=None,
pretty_name=None,
pbasename=None,
ppbasename=None,
variables=None,
**kwds,
):
if MinMaxDerivativeStatisticsBase in cls.__mro__:
return super().__new__(
cls,
name=name,
pretty_name=pretty_name,
variables=variables,
output_params=None,
**kwds,
)
else:
return super().__new__(
cls,
name=name,
pretty_name=pretty_name,
input_fields=None,
output_params=None,
**kwds,
)
@debug
def __init__(
self,
field,
components=None,
coeffs=None,
Fmin=None,
Fmax=None,
Finf=None,
all_quiet=None,
name=None,
pretty_name=None,
pbasename=None,
ppbasename=None,
variables=None,
**kwds,
):
"""
Initialize a MinMaxField Statistics operator frontend.
Available operator backends are PYTHON and OPENCL.
MinMaxFieldStatistics can compute some commonly required Field statistics:
Fmin: component-wise min values of the field.
Fmax: component-wise max values of the field.
Finf: component-wise max values of the absolute value of
the field (computed using Fmin and Fmax).
All statistics are only computed if explicitely requested by user,
unless required to compute another user-required statistic, see Notes.
All statistics may also be additionaly scaled by a coefficient.
Parameters
----------
field: Field
The continuous field on which the statistics will be computed.
components: array like of ints, optional
The components on which the statistics are computed,
defaults to all components.
coeffs: dict of array like of coefficients, optional
Optional scaling of the statistics.
Scaling factor should be a scalar or an array-like of scalars for
each components. If not given, defaults to 1 for all statistics.
F...: TensorParameter or boolean, optional
At least one statistic should be specified (either by boolean or TensorParameter).
TensorParameters should be of shape (nb_components,), see Notes.
If set to True, the TensorParameter will be generated automatically.
Autogenerated TensorParameters that are not required by the user
(ie. left to None or False during __init__) are, if required,
generated but set to be quiet.
all_quiet: bool
Set all autogenerated TensorParameter to be quiet.
name: str, optional
Name of this operator.
pretty_name: str, optional
Pretty name of this operator.
pbasename: str, optional
Parameters basename for created parameters.
Defaults to field.name.
ppbasename: str, optional
Parameters pretty basename for created parameters.
Defaults to pbasename.
variables: dict
Dictionary of fields as keys and topologies as values.
implementation: hysop.constants.Implementation, optional
base_kwds: dict, optional
Base class keyword arguments as a dictionnary.
kwds:
Extra keyword arguments passed towards operator backend implementation.
Attributes:
-----------
F...: TensorParameter
All generated tensor parameters.
Unused statistics are set to None.
Notes
-----
nb_components = min(field.nb_components, len(components)).
About statistics:
Finf requires to compute Fmin and Fmax and will have value:
Finf = Sinf * max( abs(Smin*Fmin), abs(Smax*Fmax))
where Sinf, Smin and Smax are the scaling coefficients defined in coeffs.
"""
components = to_tuple(first_not_None(components, range(field.nb_components)))
check_instance(field, Field)
check_instance(
components,
tuple,
values=int,
allow_none=True,
minval=0,
maxval=field.nb_components - 1,
)
check_instance(
coeffs, dict, keys=str, values=(int, float, np.number), allow_none=True
)
check_instance(
variables,
dict,
keys=Field,
values=CartesianTopologyDescriptors,
allow_none=True,
)
check_instance(name, str, allow_none=True)
check_instance(pbasename, str, allow_none=True)
check_instance(ppbasename, str, allow_none=True)
check_instance(all_quiet, bool, allow_none=True)
coeffs = first_not_None(coeffs, {})
name = first_not_None(name, f"MinMax[{field.name}]")
pretty_name = first_not_None(pretty_name, "|{}|∞".format(field.pretty_name))
variables = first_not_None(variables, {field: None})
all_quiet = first_not_None(all_quiet, False)
parameters = self.build_parameters(
field=field,
components=components,
all_quiet=all_quiet,
Fmin=Fmin,
Fmax=Fmax,
Finf=Finf,
pbasename=pbasename,
ppbasename=ppbasename,
)
output_params = set(p for p in parameters.values() if (p is not None))
if MinMaxDerivativeStatisticsBase in self.__class__.__mro__:
super().__init__(
name=name,
pretty_name=pretty_name,
variables=variables,
output_params=output_params,
**kwds,
)
self.has_derivative = True
else:
input_fields = {field: variables[field]}
super().__init__(
name=name,
pretty_name=pretty_name,
input_fields=input_fields,
output_params=output_params,
**kwds,
)
self.has_derivative = False
for pname, param in parameters.items():
setattr(self, pname, param)
coeffs.setdefault(pname, 1)
unused_coeffs = set(coeffs.keys()) - set(parameters.keys())
if unused_coeffs:
msg = "The following coefficients are not needed: {}"
msg = msg.format(unused_coeffs)
raise ValueError(unused_coeffs)
self._components = components
self._parameters = parameters
self._coeffs = coeffs
self._field = field
[docs]
def discretize(self):
super().discretize()
if self.has_derivative:
self._dfield = self.get_output_discrete_field(self._field)
else:
self._dfield = self.get_input_discrete_field(self._field)
for param in self._parameters.values():
if param is not None:
param.min_max_dfield = self._dfield
# Collect values from all MPI process
if self.mpi_params.size == 1:
self._collect_min = lambda e: e
self._collect_max = lambda e: e
else:
comm = self.mpi_params.comm
Fmin, Fmax = self.Fmin, self.Fmax
if self.Fmax is not None:
sendbuff = np.zeros_like(Fmax.tensor_value)
recvbuff = np.zeros_like(Fmax.tensor_value)
def _collect_max(val, sendbuff=sendbuff, recvbuff=recvbuff):
sendbuff[...] = val
comm.Allreduce(sendbuff, recvbuff, op=MPI.MAX)
return recvbuff.copy()
else:
_collect_max = None
if self.Fmin is not None:
sendbuff = np.zeros_like(Fmin.tensor_value)
recvbuff = np.zeros_like(Fmin.tensor_value)
def _collect_min(val, sendbuff=sendbuff, recvbuff=recvbuff):
sendbuff[...] = val
comm.Allreduce(sendbuff, recvbuff, op=MPI.MIN)
return recvbuff.copy()
else:
_collect_min = None
self._collect_max = _collect_max
self._collect_min = _collect_min
[docs]
def compute_statistics(self, **kwds):
"""Backend agnostic computation of min and max parameters."""
dfield, components, coeffs = self._dfield, self._components, self._coeffs
Fmin, Fmax, Finf = self.Fmin, self.Fmax, self.Finf
# For now field.min/field.max take into account ghosts...
# This is because pyopencl.reduction does not support reductions on array views.
# Here we force synchronize ghosts in every direction, including diagonals,
# in order to overwrite potential bad boundary values.
dfield.exchange_ghosts()
if Fmin is not None:
fmin = Fmin.tensor_value
for i in components:
fmin[i] = dfield.data[i].min().get()
Fmin.value = self._collect_min(fmin * coeffs["Fmin"])
if Fmax is not None:
fmax = Fmax.tensor_value
for i in components:
fmax[i] = dfield.data[i].max().get()
Fmax.value = self._collect_max(fmax * coeffs["Fmax"])
if Finf is not None:
self.Finf.value = (
np.maximum(np.abs(Fmin()), np.abs(Fmax())) * coeffs["Finf"]
)
[docs]
class MinMaxDerivativeStatisticsBase(MinMaxFieldStatisticsBase):
"""
Abstract operator base to compute min and max statistics on the derivative
of a specific scalar field.
"""
@debug
def __new__(
cls,
F,
dF=None,
A=None,
derivative=None,
direction=None,
Fmin=None,
Fmax=None,
Finf=None,
coeffs=None,
all_quiet=False,
name=None,
pretty_name=None,
pbasename=None,
ppbasename=None,
variables=None,
**kwds,
):
return super().__new__(
cls,
field=dF,
coeffs=coeffs,
Fmin=Fmin,
Fmax=Fmax,
Finf=Finf,
name=name,
pretty_name=pretty_name,
pbasename=pbasename,
variables=variables,
F=F,
dF=dF,
A=A,
derivative=derivative,
direction=direction,
**kwds,
)
@debug
def __init__(
self,
F,
dF=None,
A=None,
derivative=None,
direction=None,
Fmin=None,
Fmax=None,
Finf=None,
coeffs=None,
all_quiet=False,
name=None,
pretty_name=None,
pbasename=None,
ppbasename=None,
variables=None,
**kwds,
):
"""
Initialize an MinMaxDerivativeStatisticsBase.
MinMaxDerivativeStatistics can compute some commonly required Field derivative statistics:
Fmin: min value of a derivative of the field.
Fmax: max value of a derivative of the field.
Finf: max value of the absolute value of a
derivative of the field (computed using Fmin and Fmax).
First compute the derivative of a scalar field F in a given direction
at a given order and on a given backend out of place in a specific output scalar field dF.
The derivative is then possibly scaled by another field/parameter/value A.
After the scaled derivative has been computed, compute user requested statistics
(min and max values) on this new field and scale those statistics by other scaling
parameters stored in coeffs.
1) Compute derivative
dF = alpha * d^n(F)/dXj**n
2) Compute statistics
Fmin = Smin * min(dF)
Fmax = Smax * max(dF)
Finf = Sinf * max(|Fmin|, |Fmax|)
where F is an input field
dF is an output field (by default a temporary field).
n = derivative order > 0
alpha = A, where A is a Field, a Parameter or a scalar.
Fmin = created or supplied TensorParameter.
Fmax = created or supplied TensorParameter.
Finf = created or supplied TensorParameter.
Smin = coeffs['Fmin']
Smax = coeffs['Fmax']
Sinf = coeffs['Finf']
Statistics are only computed if explicitely requested by user,
unless required to compute another user-required statistic, see Notes.
Parameters
----------
F: hysop.field.continuous_field.Field
Continuous field as input.
dF: hysop.field.continuous_field.Field, optional
Continuous field to be written.
Some backend may allow inplace differentiation.
By default a temporary scalar field is created,
which means that the computation of the derivative
may be discarded after this operator has been applied
because of temporary buffer sharing between operators.
A: numerical value, ScalarParameter or Field, optional
Scaling field/parameter/value for convenience.
Defaults to no scaling.
derivative: int, optional
Which derivative to generate.
Defaults to 1.
direction: int, optional
Directions in which to take the derivative.
Defaults to 0.
F...: TensorParameter or boolean, optional
The output parameters that will contain the statistics.
At least one statistic should be specified (either by boolean or TensorParameter).
TensorParameters should be of shape (1,).
If set to True, the TensorParameter will be generated automatically.
Autogenerated TensorParameters that are not required by the user
(ie. left to None or False during __init__) are, if required,
generated but set to be quiet.
Autogenerated parameters can be retrieved using the 'Fmin', 'Fmax'
and 'Finf' attributes.
all_quiet: bool, optional defaults to False
Set all autogenerated TensorParameter views to be quiet.
coeffs: dict of array like of coefficients, optional
Optional scaling of the statistics.
Scaling factor should be a scalar or an array-like of scalars for
each components. If not given, defaults to 1 for all statistics.
name: str, optional
Name of this operator.
pretty_name: str, optional
Pretty name of this operator.
pbasename: str, optional
Parameters basename for created parameters.
Defaults to field.name.
ppbasename: str, optional
Parameters pretty basename for created parameters.
Defaults to pbasename.
variables: dict, optional
Dictionary of fields as keys and topologies as values.
implementation: hysop.constants.Implementation, optional
Specify underlying backend implementation.
Target implementation, should be contained in available_implementations().
If None, implementation will be set to default_implementation().
base_kwds: dict, optional
Base class keyword arguments as a dictionnary.
kwds:
Extra keyword arguments passed towards operator backend implementation.
Attributes:
-----------
F...: TensorParameter
All generated tensor parameters.
Unused statistics are set to None.
Notes
-----
About statistics:
Finf requires to compute Fmin and Fmax and will have value:
Finf = Sinf * max( abs(Smin*Fmin), abs(Smax*Fmax))
where Sinf, Smin and Smax are the scaling coefficients defined in coeffs.
"""
if variables is None:
msg = "The variables parameter cannot be None."
raise ValueError(msg)
assert F in variables, variables
if dF is None:
dF = F.tmp_like("dF", nb_components=1)
variables.setdefault(dF, variables[F])
super().__init__(
field=dF,
coeffs=coeffs,
Fmin=Fmin,
Fmax=Fmax,
Finf=Finf,
name=name,
pretty_name=pretty_name,
pbasename=pbasename,
variables=variables,
F=F,
dF=dF,
A=A,
derivative=derivative,
direction=direction,
**kwds,
)